Source Visibility with agilepy¶

The AGILE satellite operated from April 23, 2007 to January 18, 2024 observing thousands of astrophysical sources.

The AGILE spacecraft operated in "pointing mode" from the beginning of the mission to October 15, 2009, completing 101 observation blocks (OBs). The OBs usually consisted of predefined long exposures, drifting about 1 deg per day with respect to the initial boresight direction to obey solar panels constraints.

In November 2009, the attitude control system was reconfigured, and scientific operations were performed "spinning mode" until the end of the mission. AGILE scanned ≈ 80% sky daily (exposure of ≈ 7 · 106 cm2 s) with an angular velocity of about 0.8 deg s−1, performing 200 passes per day on the same sky region.

The visibility of a source thus depend on the off-axis angle with respect to the satellite pointing direction.

The purpose of this notebook is to show how to use the class AGVisibility to determine the AGILE pointing direction and the offset of a source of interest. The class requires only the "log files", which are the AGILE spacecraft files, to perform its task. If a Fermi Spacecraft file is provided, the same operation can be performed on the file for comparison purposes.

In [1]:
# Import the relevant class
from agilepy.api.AGVisibility import AGVisibility
from agilepy.utils.AstroUtils import AstroUtils

from astropy.coordinates import SkyCoord
from astropy.time import Time

Setup¶

This tutorial uses the available test dataset on blazar 3C 454.3, from 2010-11-13 to 2010-11-21 (MJD 55513-55521).

In [2]:
# Prepare the YAML Configuration File.
# You can create it manually or use AGVisibility.getConfiguration()

outputDir    = "/home/flareadvocate/workspace/shared_dir/"
confFilePath = "/home/flareadvocate/workspace/shared_dir/tutorial_visibility.yaml"
logFile      = "$AGILE/agilepy-test-data/test_dataset_agn/LOG/LOG.index"
fermi_SCFile = "$AGILE/agilepy-test-data/visibility/Fermi_test_SC00.fits"

AGVisibility.getConfiguration(
    confFilePath=confFilePath,
    # Output
    outputDir=outputDir,
    userName="my_name",
    sourceName="vis-source",
    verboselvl=0,
    # Input
    logFile=logFile,
    fermiSpacecraftFile=fermi_SCFile,
    # Selection
    step=1,
    timeType="tt",
    tmin="null",
    tmax="null",
    # Source
    coord1="null",
    coord2="null",
    frame="icrs",
)

The YAML file has the following sections and field:

  • output: common to all agilepy classes, sets the following output fields: outdir, filenameprefix, sourcename, username, verboselvl
  • input: input spacecraft files.
    • logFile: path to the index file listing relevant AGILE log files to be extracted.
    • fermiSpacecraftFile: the path to the Fermi spacecraft file. Optional.
  • selection: arguments that define the timing of the information extracted.
    • step: time interval in seconds between two consecutive points in the resulting table. Minimum accepted value: 0.1.
    • timeType: the time format of tmin and tmax. AGILE seconds can be expressed as tt.
    • tmin: inferior observation time limit to analize. If not set, default value is TT = 104371200 (2007-04-23 00:00:00)
    • tmax: superior observation time limit to analize. If not set, default value is TT = 632620800 (2024-01-18 00:00:00)
  • source: optional, include arguments to define the coordinates of the source used to compute the offset. The coordinates can also be set or changed when using the class. They are defined by:
    • frame: reference frame for astropy.coordinates.SkyCoord
    • coord1, coord2: source coordinates in degrees.
In [3]:
# Create the AGVisibility object
ag_vis = AGVisibility(confFilePath)

# Print Options
ag_vis.printOptions()
Log level set to WARNING and output to /home/flareadvocate/workspace/shared_dir/my_name_vis-source_20251005-231634/logs
{ 'input': { 'fermi_spacecraft_file': '/home/flareadvocate/agiletools/agilepy-test-data/visibility/Fermi_test_SC00.fits',
             'logfile': '/home/flareadvocate/agiletools/agilepy-test-data/test_dataset_agn/LOG/LOG.index'},
  'output': { 'filenameprefix': 'visibility_product',
              'outdir': PosixPath('/home/flareadvocate/workspace/shared_dir/my_name_vis-source_20251005-231634'),
              'sourcename': 'vis-source',
              'username': 'my_name',
              'verboselvl': 0},
  'selection': { 'step': 1,
                 'timetype': 'TT',
                 'tmax': 632620800,
                 'tmin': 104371200},
  'source': {'coord1': None, 'coord2': None, 'frame': 'icrs'}}
In [4]:
# You can use the setOptions() method to change the configuration:

# Set coordinates of 3C 454.3
ag_vis.setOptions(coord1=343.496, coord2=16.151, frame='icrs')

# Set time range appropriate for the test dataset
ag_vis.setOptions(tmin=AstroUtils.convert_time_to_agile_seconds(Time(55513, format='mjd')),
                  tmax=AstroUtils.convert_time_to_agile_seconds(Time(55521, format='mjd')),
                  timetype="TT"
                  )

# Check configuration
ag_vis.printOptions()
{ 'input': { 'fermi_spacecraft_file': '/home/flareadvocate/agiletools/agilepy-test-data/visibility/Fermi_test_SC00.fits',
             'logfile': '/home/flareadvocate/agiletools/agilepy-test-data/test_dataset_agn/LOG/LOG.index'},
  'output': { 'filenameprefix': 'visibility_product',
              'outdir': PosixPath('/home/flareadvocate/workspace/shared_dir/my_name_vis-source_20251005-231634'),
              'sourcename': 'vis-source',
              'username': 'my_name',
              'verboselvl': 0},
  'selection': { 'step': 1,
                 'timetype': 'TT',
                 'tmax': 217382400.0,
                 'tmin': 216691200.0},
  'source': {'coord1': 343.496, 'coord2': 16.151, 'frame': 'icrs'}}

Note on time conversion: the agilepy.utils.AstroUtils.AstroUtils class offers utility functions for the conversion between AGILE seconds and other data formats, handled through the astropy.time.Time class.

In [5]:
# Convert from AGILE seconds to astropy Time
t_example = AstroUtils.convert_time_from_agile_seconds(104371200)
display(f"TT=104371200 corresponds to MJD {t_example.mjd}, ISO {t_example.iso}, UNIX {t_example.unix}")

# Convert from astropy Time to AGILE seconds
t_example = AstroUtils.convert_time_to_agile_seconds(Time(54213, format='mjd'))
display(f"MJD=54213.0 corresponds to TT {t_example}")
'TT=104371200 corresponds to MJD 54213.0, ISO 2007-04-23 00:00:00.000, UNIX 1177286400.0'
'MJD=54213.0 corresponds to TT 104371200.0'

AGILE Visibility¶

The AGILE pointing coordinates can be computed with AGVisibility.computePointingDirection(). The arguments of the function are taken by default from the configuration, or they can be overriden if provided explicitly as arguments of the method.

The AGILE data is extracted from the columns TIME, ATTITUDE_RA_Y, ATTITUDE_DEC_Y of the log files and stored in an astropy.table.Table.

The output columns provide the start and stop times of the bin analysed in TT and MJD, the AGILE pointing coordinates (Right Ascension and Declination) and the mean offset of the source in degrees in the bin. The distance between two consecutive bins is the step parameter.

The table is saved by default in the agilepy output directory as a .csv file.

In [6]:
agile_visibility_table = ag_vis.computePointingDirection(writeFiles=True)

# Show
display(agile_visibility_table)

# Note: The table can also be accessed as a property of AGVisibility:
# ag_vis.agileVisibility
Table length=691052
TT_startTT_stopAGILE_RAAGILE_DECMJD_startMJD_stopoffaxis_angle_deg
float64float64float64float64float64float64float64
216691200.0216691200.1222.2467651367187572.388793945312555513.055513.0000011574183.43257364718369
216691201.0216691201.1219.7959136962890672.3218383789062555513.0000115740855513.0000127314884.07688376120517
216691202.0216691202.1217.3642425537109472.2247390747070355513.00002314814655513.00002430555584.72213619894937
216691203.0216691203.1214.9628295898437572.0969085693359455513.0000347222255513.0000358796385.36813315305827
216691204.0216691204.1212.596603393554771.9396362304687555513.000046296355513.000047453786.01495135877619
216691205.0216691205.1210.2767486572265671.7533798217773455513.0000578703755513.0000590277886.66169258017244
216691206.0216691206.1208.003952026367271.5381546020507855513.00006944444655513.0000706018587.31008697219828
216691207.0216691207.1205.7929840087890671.2930984497070355513.00008101851655513.00008217592487.95946833316754
.....................
217382392.1217382392.2198.8926086425781265.1190261840820355520.99990856481655520.99990972222494.42071818104152
217382393.1217382393.2197.4219360351562564.6984939575195355520.99992013888555520.99992129629495.11427439693121
217382394.1217382394.2195.9965972900390664.2629089355468855520.9999317129655520.9999328703795.80856711789367
217382395.1217382395.2194.616500854492263.8131523132324255520.9999432870455520.9999444444596.50327748643129
217382396.1217382396.2193.2785034179687563.3522796630859455520.9999548611155520.9999560185297.1967921444698
217382397.1217382397.2191.9813385009765662.87891006469726655520.99996643518555520.99996759259497.89085242403705
217382398.1217382398.2190.7301940917968862.3926277160644555520.99997800925555520.9999791666798.58504254443122
217382399.1217382399.2189.5220336914062561.894519805908255520.9999895833355520.9999907407499.27936207655596

Fermi Visibility¶

The Fermi pointing coordinates can be computed with AGVisibility.getFermiPointing(). The arguments of the function are taken by default from the configuration, or they can be overriden if provided explicitly as arguments of the method.

The Fermi data is extracted from the columns START, STOP, RA_SCZ, DEC_SCZ of the spacecraft file and stored in an astropy.table.Table.

The output columns provide the start and stop times of the bin analysed in TT, MJD, Fermi MET seconds, the Fermi pointing coordinates (Right Ascension and Declination) and the mean offset of the source in degrees in the bin.

The table is saved by default in the agilepy output directory as a .csv file.

In [7]:
fermi_visibility_table = ag_vis.getFermiPointing(writeFiles=True)

# Show the table, accessed as a property of AGVisibility:
display(ag_vis.fermiVisibility)
WARNING: UnitsWarning: 'Gauss' did not parse as fits unit: At col 0, Unit 'Gauss' not supported by the FITS standard.  If this is meant to be a custom unit, define it with 'u.def_unit'. To have it recognized inside a file reader or other code, enable it with 'u.add_enabled_units'. For details, see https://docs.astropy.org/en/latest/units/combining_and_defining.html [astropy.units.core]
WARNING: UnitsWarning: 'Earth_Radii' did not parse as fits unit: At col 0, Unit 'Earth_Radii' not supported by the FITS standard.  If this is meant to be a custom unit, define it with 'u.def_unit'. To have it recognized inside a file reader or other code, enable it with 'u.add_enabled_units'. For details, see https://docs.astropy.org/en/latest/units/combining_and_defining.html [astropy.units.core]
Table length=19400
MET_STARTMET_STOPRA_SCZDEC_SCZTT_startTT_stopMJD_startMJD_stopoffaxis_angle_deg
ssdegdeg
float64float64float32float32float64float64float64float64float32
311299218.6311299248.6266.6657470.37467216691218.60000002216691248.6000000255513.00021527777455513.000562570.39624
311299248.6311299278.6269.5991269.69561216691248.60000002216691278.6000000255513.000562555513.0009097222269.30859
311299278.6311299308.6272.395468.99575216691278.60000002216691308.6000000255513.0009097222255513.00125694444468.210075
311299308.6311299338.6275.0628468.274506216691308.60000002216691338.6000000255513.00125694444455513.0016041666767.1011
311299338.6311299368.6277.6157267.53844216691338.60000002216691368.6000000255513.0016041666755513.0019513888965.98142
311299368.6311299398.6280.0590866.78771216691368.60000002216691398.6000000255513.0019513888955513.0022986111164.85262
311299398.6311299428.6282.403766.0241216691398.60000002216691428.6000000255513.0022986111155513.0026458333363.714382
...........................
311988676.6311988706.6321.1410225.867832217380676.60000002217380706.6000000255520.9800532407455520.9804004629622.9716
311988706.6311988736.6322.4569425.6432217380706.60000002217380736.6000000255520.9804004629655520.98074768518621.785173
311988736.6311988766.6323.7742625.435379217380736.60000002217380766.6000000255520.98074768518655520.9810949074120.606659
311988766.6311988796.6325.093525.246605217380766.60000002217380796.6000000255520.9810949074155520.98144212962519.43865
311988796.6311988824.085537326.412525.07553217380796.60000002217380824.0855370255520.98144212962555520.9817602492718.285112
311990302.64574337311990332.686.6461-48.26763217382302.64574337217382332.6000000255520.9988732146255520.999219907404110.67387
311990332.6311990362.688.04295-47.49842217382332.60000002217382362.6000000255520.99921990740455520.99956712963111.59731
311990362.6311990392.689.429375-46.734695217382362.60000002217382392.6000000255520.9995671296355520.99991435185112.53722

Plot¶

AGVisibility offers an utility function to plot the time series of the Offaxis angle. Parameters are:

  • sourceCoordinates (astropy.coordinates.SkyCoord): the coordinates of the source used to compute the offaxis angle, If None, the offaxis angle must already be computed in the visibility table.
  • maxOffaxis (float): maximum offaxis angle in degrees to be shown in the visibility plot.
  • mjd_limits (tuple(float,float), optional): Plot limits in MJD.
  • plotFermi (bool): if True, plot the Fermi offaxis angle.
  • plotHistogram (bool): if True, plot the histogram of the offaxis angles.
  • saveImage (bool): if True, save the images in the output directory.
  • showPositionPanels (bool): if True, show the RA and DEC time series of the satellite direction.
In [8]:
# Plot
ag_vis.plotVisibility(mjd_limits=(55519, 55519.5), maxOffaxis=60,
                      plotFermi=True, plotHistogram=True, saveImages=True,
                      showPositionPanels=False
                      )
Out[8]:
[PosixPath('/home/flareadvocate/workspace/shared_dir/my_name_vis-source_20251005-231634/plots/offaxis_ra343.4960_+16.1510.reg_MJD_55519_55519.5.png'),
 PosixPath('/home/flareadvocate/workspace/shared_dir/my_name_vis-source_20251005-231634/plots/offaxis_hist_ra343.4960_+16.1510.reg_MJD_55519_55519.5.png')]
In [ ]: